# # Contracting cube
#
# In this demo we simulate a unit cube that is fixed
# at $x = 0$ and free at $x = 1$. We use
# a transversally isotropic material with fiber oriented
# in the $x$-direction.
import dolfin
try:
from dolfin_adjoint import (
Constant,
DirichletBC,
Expression,
UnitCubeMesh,
interpolate,
Mesh,
)
except ImportError:
from dolfin import (
Constant,
DirichletBC,
interpolate,
Expression,
UnitCubeMesh,
Mesh,
)
import pulse
from fenics_plotly import plot
pulse.iterate.logger.setLevel(10)
# Create mesh
N = 6
mesh = UnitCubeMesh(N, N, N)
# Create subdomains
class Free(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] > (1.0 - dolfin.DOLFIN_EPS) and on_boundary
class Fixed(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] < dolfin.DOLFIN_EPS and on_boundary
# Create a facet fuction in order to mark the subdomains
ffun = dolfin.MeshFunction("size_t", mesh, 2)
ffun.set_all(0)
# Mark the first subdomain with value 1
fixed = Fixed()
fixed_marker = 1
fixed.mark(ffun, fixed_marker)
# Mark the second subdomain with value 2
free = Free()
free_marker = 2
free.mark(ffun, free_marker)
# Create a cell function (but we are not using it)
cfun = dolfin.MeshFunction("size_t", mesh, 3)
cfun.set_all(0)
# Collect the functions containing the markers
marker_functions = pulse.MarkerFunctions(ffun=ffun, cfun=cfun)
# Create mictrotructure
V_f = dolfin.VectorFunctionSpace(mesh, "CG", 1)
# Fibers
f0 = interpolate(Expression(("1.0", "0.0", "0.0"), degree=1), V_f)
# Sheets
s0 = interpolate(Expression(("0.0", "1.0", "0.0"), degree=1), V_f)
# Fiber-sheet normal
n0 = interpolate(Expression(("0.0", "0.0", "1.0"), degree=1), V_f)
# Collect the mictrotructure
microstructure = pulse.Microstructure(f0=f0, s0=s0, n0=n0)
# Create the geometry
geometry = pulse.Geometry(
mesh=mesh,
marker_functions=marker_functions,
microstructure=microstructure,
)
# Use the default material parameters
material_parameters = pulse.HolzapfelOgden.default_parameters()
# Select model for active contraction
active_model = pulse.ActiveModels.active_strain
# active_model = "active_stress"
# Set the activation
activation = Constant(0.0)
# Create material
material = pulse.HolzapfelOgden(
active_model=active_model,
parameters=material_parameters,
activation=activation,
)
# Make Dirichlet boundary conditions
def dirichlet_bc(W):
V = W if W.sub(0).num_sub_spaces() == 0 else W.sub(0)
return DirichletBC(V, Constant((0.0, 0.0, 0.0)), fixed)
# Make Neumann boundary conditions
neumann_bc = pulse.NeumannBC(traction=Constant(0.0), marker=free_marker)
# Collect Boundary Conditions
bcs = pulse.BoundaryConditions(dirichlet=(dirichlet_bc,), neumann=(neumann_bc,))
# Create problem
problem = pulse.MechanicsProblem(geometry, material, bcs)
# Solve problem
pulse.iterate.iterate(problem, activation, 0.1)
2021-11-23 07:14:57,112 [630] DEBUG pulse.iterate: Control: [0.0]
2021-11-23 07:14:57,115 [630] DEBUG pulse.iterate: Target: [0.1]
2021-11-23 07:14:57,117 [630] DEBUG pulse.iterate: Intial number of steps: 5 with step size 0.02
2021-11-23 07:14:57,118 [630] INFO pulse.iterate: Iterating to target control (f_21)...
2021-11-23 07:14:57,119 [630] INFO pulse.iterate: Current control: f_21 = 0.000
2021-11-23 07:14:57,120 [630] INFO pulse.iterate: Target: 0.100
2021-11-23 07:14:57,120 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:14:57,121 [630] DEBUG pulse.iterate: Maximum difference: 1.000e-01
2021-11-23 07:14:57,122 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:14:57,123 [630] DEBUG pulse.iterate: Current control: 0.020
2021-11-23 07:14:59,116 [630] DEBUG pulse.iterate: Solver did not converge...
2021-11-23 07:14:59,117 [630] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-23 07:14:59,117 [630] DEBUG pulse.iterate: Reduce control step
2021-11-23 07:14:59,118 [630] DEBUG pulse.iterate: Assign old state
2021-11-23 07:14:59,124 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:14:59,125 [630] DEBUG pulse.iterate: Maximum difference: 1.000e-01
2021-11-23 07:14:59,126 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:14:59,126 [630] DEBUG pulse.iterate: Current control: 0.010
2021-11-23 07:15:05,130 [630] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 07:15:05,131 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:05,131 [630] DEBUG pulse.iterate: Maximum difference: 9.000e-02
2021-11-23 07:15:05,133 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:05,134 [630] DEBUG pulse.iterate: Maximum difference: 9.000e-02
2021-11-23 07:15:05,134 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:05,135 [630] DEBUG pulse.iterate: Current control: 0.020
2021-11-23 07:15:07,608 [630] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 07:15:07,609 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:07,610 [630] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2021-11-23 07:15:07,611 [630] DEBUG pulse.iterate: Adapt step size. New step size: 0.015
2021-11-23 07:15:07,612 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:07,614 [630] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2021-11-23 07:15:07,615 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:07,615 [630] DEBUG pulse.iterate: Current control: 0.035
2021-11-23 07:15:11,144 [630] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 07:15:11,145 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:11,146 [630] DEBUG pulse.iterate: Maximum difference: 6.500e-02
2021-11-23 07:15:11,147 [630] DEBUG pulse.iterate: Adapt step size. New step size: 0.022
2021-11-23 07:15:11,148 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:11,149 [630] DEBUG pulse.iterate: Maximum difference: 6.500e-02
2021-11-23 07:15:11,150 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:11,151 [630] DEBUG pulse.iterate: Current control: 0.058
2021-11-23 07:15:14,163 [630] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 07:15:14,164 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:14,165 [630] DEBUG pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 07:15:14,166 [630] DEBUG pulse.iterate: Adapt step size. New step size: 0.034
2021-11-23 07:15:14,167 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:14,168 [630] DEBUG pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 07:15:14,169 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:14,170 [630] DEBUG pulse.iterate: Current control: 0.091
2021-11-23 07:15:15,794 [630] DEBUG pulse.iterate: Solver did not converge...
2021-11-23 07:15:15,795 [630] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-23 07:15:15,795 [630] DEBUG pulse.iterate: Reduce control step
2021-11-23 07:15:15,796 [630] DEBUG pulse.iterate: Assign old state
2021-11-23 07:15:15,803 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:15,804 [630] DEBUG pulse.iterate: Maximum difference: 4.250e-02
2021-11-23 07:15:15,806 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:15,807 [630] DEBUG pulse.iterate: Current control: 0.074
2021-11-23 07:15:18,820 [630] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 07:15:18,822 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:18,822 [630] DEBUG pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 07:15:18,823 [630] DEBUG pulse.iterate: Adapt step size. New step size: 0.025
2021-11-23 07:15:18,824 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:18,826 [630] DEBUG pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 07:15:18,827 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:18,828 [630] DEBUG pulse.iterate: Current control: 0.100
2021-11-23 07:15:20,409 [630] DEBUG pulse.iterate: Solver did not converge...
2021-11-23 07:15:20,411 [630] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-23 07:15:20,411 [630] DEBUG pulse.iterate: Reduce control step
2021-11-23 07:15:20,412 [630] DEBUG pulse.iterate: Assign old state
2021-11-23 07:15:20,419 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:20,420 [630] DEBUG pulse.iterate: Maximum difference: 2.563e-02
2021-11-23 07:15:20,421 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:20,422 [630] DEBUG pulse.iterate: Current control: 0.087
2021-11-23 07:15:22,909 [630] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 07:15:22,911 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:22,912 [630] DEBUG pulse.iterate: Maximum difference: 1.297e-02
2021-11-23 07:15:22,913 [630] DEBUG pulse.iterate: Adapt step size. New step size: 0.019
2021-11-23 07:15:22,914 [630] DEBUG pulse.iterate: Check target reached: NO
2021-11-23 07:15:22,915 [630] DEBUG pulse.iterate: Maximum difference: 1.297e-02
2021-11-23 07:15:22,916 [630] DEBUG pulse.iterate: Change step size for final iteration
2021-11-23 07:15:22,916 [630] DEBUG pulse.iterate: Try new control
2021-11-23 07:15:22,917 [630] DEBUG pulse.iterate: Current control: 0.100
2021-11-23 07:15:24,954 [630] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-23 07:15:24,955 [630] DEBUG pulse.iterate: Check target reached: YES!
2021-11-23 07:15:24,956 [630] DEBUG pulse.iterate: Check target reached: YES!
([Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 47),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 195),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 246),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 313),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 372),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 466),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 552),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 595)],
[Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 45),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 193),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 244),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 311),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 370),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 464),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 550),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 593)])
# Get displacement and hydrostatic pressure
u, p = problem.state.split(deepcopy=True)
V = dolfin.VectorFunctionSpace(mesh, "CG", 1)
u_int = interpolate(u, V)
new_mesh = Mesh(mesh)
dolfin.ALE.move(new_mesh, u_int)
fig = plot(mesh, show=False)
fig.add_plot(plot(new_mesh, color="red", show=False))
fig.show()